Introduction

Pymatgen supports reading of most common file formats, including the Crystallographic Information File and various input and output files of computational codes like VASP. However, it is often easier and quicker to directly query for structures from online sources. Though private databases such as the Inorganic Crystal Structure Database are not open, there are open sources such as the Materials Project and the Crystallographic Open Database (COD) where one can obtain crystal structures.

Written using:

  • pymatgen==2018.3.13

Materials Project

Pymatgen contains a high-level interface to the Materials Project, which can be used to query for structures very easily.


In [1]:
from pymatgen.ext.matproj import MPRester
# Note that you will need to add your Materials API key in your .pmgrc.yaml file as "PMG_MAPI_KEY". 
# Alternatively, you will need to supply the API key as an arg to MPRester.
mpr = MPRester()

In [2]:
# Querying by formula only.
structures = mpr.get_structures("Li2O")
print(structures)


[Structure Summary
Lattice
    abc : 3.2910717923597561 3.2910718996250861 3.2910720568557887
 angles : 60.129710432884849 60.129709521376753 60.129703130390972
 volume : 25.279668381289056
      A : 2.91738857 0.097894369999999994 1.5200046599999999
      B : 0.96463405999999996 2.7550356100000002 1.5200046599999999
      C : 0.13320635 0.097894430000000004 3.28691771
PeriodicSite: O (0.0000, 0.0000, 0.0000) [0.0000, 0.0000, 0.0000]
PeriodicSite: Li (3.0121, 2.2136, 4.7463) [0.7502, 0.7502, 0.7502]
PeriodicSite: Li (1.0031, 0.7372, 1.5806) [0.2498, 0.2498, 0.2498], Structure Summary
Lattice
    abc : 5.1517948200000001 3.1404278300000001 5.9334081599999999
 angles : 90.0 90.0 90.0
 volume : 95.995660249910003
      A : 5.1517948200000001 0.0 0.0
      B : 0.0 3.1404278300000001 0.0
      C : 0.0 0.0 5.9334081599999999
PeriodicSite: Li (0.0631, 0.7851, 0.9438) [0.0122, 0.2500, 0.1591]
PeriodicSite: Li (0.7268, 0.7851, 3.4367) [0.1411, 0.2500, 0.5792]
PeriodicSite: Li (1.8490, 2.3553, 0.4700) [0.3589, 0.7500, 0.0792]
PeriodicSite: Li (2.5128, 2.3553, 3.9105) [0.4878, 0.7500, 0.6591]
PeriodicSite: Li (2.6390, 0.7851, 2.0229) [0.5122, 0.2500, 0.3409]
PeriodicSite: Li (3.3027, 0.7851, 5.4634) [0.6411, 0.2500, 0.9208]
PeriodicSite: Li (4.4249, 2.3553, 2.4967) [0.8589, 0.7500, 0.4208]
PeriodicSite: Li (5.0887, 2.3553, 4.9896) [0.9878, 0.7500, 0.8409]
PeriodicSite: O (1.2677, 2.3553, 2.3664) [0.2461, 0.7500, 0.3988]
PeriodicSite: O (1.3082, 0.7851, 5.3331) [0.2539, 0.2500, 0.8988]
PeriodicSite: O (3.8436, 2.3553, 0.6003) [0.7461, 0.7500, 0.1012]
PeriodicSite: O (3.8841, 0.7851, 3.5670) [0.7539, 0.2500, 0.6012]]

In [3]:
# Querying by chemical system only.
structures = mpr.get_structures("Li-O")
for s in structures:
    print(s.formula)
# A number of Li-O structures are returned with different Li:O ratios.


Li2 O4
Li2 O1
Li8 O4
Li4 O4
Li1 O3
Li16 O16

COD

To obtain structures from COD by the COD id, you just need to know the id. However, most sophisticated searches require that you have installed mysql given that the COD does not support a proper REST API at this time.


In [4]:
from pymatgen.ext.cod import COD
cod = COD()

In [5]:
structure = cod.get_structure_by_id(1010064)
print(structure)


Full Formula (Li8 O4)
Reduced Formula: Li2O
abc   :   4.610000   4.610000   4.610000
angles:  90.000000  90.000000  90.000000
Sites (12)
  #  SP       a     b     c
---  ----  ----  ----  ----
  0  Li+   0.25  0.25  0.25
  1  Li+   0.25  0.75  0.75
  2  Li+   0.75  0.25  0.75
  3  Li+   0.75  0.75  0.25
  4  Li+   0.75  0.75  0.75
  5  Li+   0.75  0.25  0.25
  6  Li+   0.25  0.75  0.25
  7  Li+   0.25  0.25  0.75
  8  O2-   0     0     0
  9  O2-   0     0.5   0.5
 10  O2-   0.5   0     0.5
 11  O2-   0.5   0.5   0

In [6]:
structures = cod.get_structure_by_formula("Li2O")
for d in structures:
    print("COD ID: %d, Formula: %s, Spacegroup: %s" % (d["cod_id"], d["structure"].formula, d["sg"]))


/Users/shyue/repos/pymatgen/pymatgen/io/cif.py:801: UserWarning: LI parsed as L
  warnings.warn("{} parsed as {}".format(sym, v))
COD ID: 1010064, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1011372, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 9009059, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 4311895, Formula: L8 O4, Spacegroup: F m -3 m
COD ID: 1514086, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1514087, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1514088, Formula: Li6 O3, Spacegroup: R -3 m :H
COD ID: 4121514, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 4121515, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1514092, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1514093, Formula: Li7.84 O4, Spacegroup: F m -3 m
COD ID: 1514094, Formula: Li7.9984 O4, Spacegroup: F m -3 m
COD ID: 1514095, Formula: Li8.0008 O4, Spacegroup: F m -3 m
COD ID: 1514096, Formula: Li8 O4, Spacegroup: F m -3 m
COD ID: 1514097, Formula: Li8.0008 O4, Spacegroup: F m -3 m
COD ID: 1514098, Formula: Li8 O4, Spacegroup: F m -3 m

In [7]:
print(structures[0]["structure"])


Full Formula (Li8 O4)
Reduced Formula: Li2O
abc   :   4.610000   4.610000   4.610000
angles:  90.000000  90.000000  90.000000
Sites (12)
  #  SP       a     b     c
---  ----  ----  ----  ----
  0  Li+   0.25  0.25  0.25
  1  Li+   0.25  0.75  0.75
  2  Li+   0.75  0.25  0.75
  3  Li+   0.75  0.75  0.25
  4  Li+   0.75  0.75  0.75
  5  Li+   0.75  0.25  0.25
  6  Li+   0.25  0.75  0.25
  7  Li+   0.25  0.25  0.75
  8  O2-   0     0     0
  9  O2-   0     0.5   0.5
 10  O2-   0.5   0     0.5
 11  O2-   0.5   0.5   0

In [ ]: